BigDFT.Systems module
This module contains the System level class of PyBigDFT. Systems are named collections of fragments, and represent a complete system for simulation.
- class System(*args, **kwargs)[source]
A system is defined as a named collection of fragments. You can manipulate a system as if it were a standard python dictionary, however it also has helper routines for performing operations on the full system.
- property centroid
Center of mass of the system
- property central_fragment
Returns the fragment whose center of mass is closest to the centroid
- Returns:
the name of the fragment. (Fragment): the fragment object
- Return type:
(str)
- get_external_potential(units='bohr', charge_offset=False)[source]
Transform the system information into a dictionary ready to be put as an external potential.
- get_k_nearest_fragments(target, k, cutoff=None, return_type='list')[source]
Given a fragment id in a system, this computes the nearest fragment.
- Parameters:
- Returns:
- the ids of the nearest fragments, or their distances
as values in the case in which ‘dict’ is provided in the return_type argument
- Return type:
- get_nearest_fragment(target)[source]
Given a fragment id in a system, this computes the nearest fragment.
- get_net_force()[source]
Returns the net force on a system in Ha/Bohr.
- Returns:
Three values which describe the net force.
- Return type:
(list)
- get_atoms(order=None)[source]
Generator enabling to iterate on system’s atoms.
This function may be useful to control the order of the atoms positions.
- Parameters:
order (list) – order of the fragments. Useful on versions of python where the systems’ dict are not ordered.
- Yields:
BigDFT.Atoms.Atom – the next atom of the system.
- get_posinp(units='angstroem', order=None)[source]
Provide the dictionary which has to be passed to the
posinp
value of therun()
method of theSystemCalculator
class instance.
- serialize(order=None, units='bohr')[source]
Transform the system in a list that can be employed for the construction of dataframes or pandas series. :param order: list of fragments to serialize in order :type order: list :param units: the units for the positions :type units: str
- Returns:
- the serialized system as well as a lookup dictionary
that contains the order of the fragment atoms in the series
- Return type:
- subsystem(fragments)[source]
Extract a subsystem.
Create a subsystem from a collection of fragments.
- to_dataframe(**kwargs)[source]
Convert the system into a dataframe, from the py:meth:~System.serialize method.
- Parameters:
**kwargs – arguments to be passed to py:meth:System.serialize method
- dataframe_slicing(df=None)[source]
Define a dictionaries of two-element lists indicating the slicing (start and end points) of the corresponding fragment in the order provided by the dataframe
- Parameters:
df (Dataframe) – associated to the system
- Returns:
dictionary of the slicings of the fragments in the dataframe
- Return type:
- to_file(filename, **kwargs)[source]
Dump the System instance into a file.
Write the system information in the file system, according to the file extension
- Parameters:
filename (str) – path of the filename to write into. Automatically determine the format according to the file extension
**kwargs – further kwyword arguments to be passed to the writing routine
- property PointParticles
~PointParticles.PointParticles` object
- Type:
Transform the system into a `py
- Type:
class
- property electrostatic_interactions
Dictionary of the Electrostatic interactions between fragments.
- property q0
Provides the global monopole of the system given as a sum of the monopoles of the atoms.
- property qcharge
The total qcharge of a system.
- rename_fragments(fragment_mapping=None)[source]
This procedure automatically names the fragments in a system.
- Parameters:
fragment_mapping (dict) – Dictionary containing list of fragments that are additionally added to each of the original system’s fragments.
- Returns:
the same system, with the automatic naming scheme.
- Return type:
- set_atom_multipoles(logfile, correct_charge=True)[source]
After a run is completed, we have a set of multipoles defined on each atom. This routine will set those values on to each atom in the system.
- Parameters:
logfile (Logfiles.Logfile) – logfile with the multipole values.
correct_charge (bool) – currently there is an inconsistency in terms of gross charge, and this corrects it.
- set_atom_forces(logfile)[source]
After a run is completed, we have the forces on each atom in the logfile. This routine will set those values to each atom in this sytem.
- Parameters:
logfile (Logfiles.Logfile) – logfile with the forces.
- compute_matching(atlist, check_matching=True)[source]
Frequently we are passed a list of atom like objects from which we need to extract data and assign it to a system. However, a system can potentially store those atoms in any order, and may not have the same set of atoms. This helper routine creates a mapping between this list view, to the dictionary view of the system class.
- Parameters:
- Returns:
- a mapping from a system to indices in the atom list. If
an atom is not in the list, an index value of -1 is assigned.
- Return type:
(dict)
- display(colordict=None, field_vals=None, cartoon=False)[source]
Display the system using the inline visualizer of Py3DMol
- examine(axs=None, view=None)[source]
Provide quick overview of the system’s quality.
This routine displays information about the coordination numbers and the bond of the system. Useful to understand if some positions of the system need optimization.
- Parameters:
axs – list of matplotlib axis in which to plot the information. If the forces are present it should be of size 3.
view (dict) – a Fragment view of the system in which to represent the fragment forces. If absent, the fragment forces are not represented.
- Returns:
Dictionary of the information about the system
- Return type:
- set_electrons_from_log(log)[source]
BigDFT uses pseudopotentials, so this will extract from a logfile the number of actual electrons modelled for each atom in the system.
- Parameters:
log (Logfiles.Logfile) – A BigDFT run logfile.
- set_logfile_info(log)[source]
Include the information of the logfile in the fragment quantities.
- Parameters:
log (Logfiles.Logfile) – A BigDFT run logfile.
- ase_potential_energy(ase_calculator)[source]
Given a ASE calculator, calculates the potential energy of the system.
- Parameters:
ase_calculator (ase.calculators.calculator.Calculator) – ASE calculator.
- Returns:
the potential energy, in Hartree
- Return type:
- reform_superunits(mapping)[source]
Creates a new system from the provided mapping
- Parameters:
mapping (dict) – dictionary of the form {newfrag: [frag1,frag2]} defining the mapping between the old fragments and the new
- Returns:
a new system from the remapping
- Return type:
- update_positions_from_dict(posinp)[source]
Update the atomic positions of a system from a posinp dictionary.
This method only works if the order of atoms match.
- Parameters:
posinp (dict) – a posinp dictionary.
- adjacency_matrix()[source]
Generates a sparse adjacency matrix based on the connectivity of a system.
- Returns:
the sparse adjacency matrix.
- Return type:
(scipy.sparse.dok_matrix)
- fragment_view(purities, bond_orders, view=None)[source]
Returns the fragment view of the system, according to a mapping.
- Parameters:
view (dict) – Mapping of the fragments to be taken initially.
btool (BigDFT.PostProcessing.BigDFTool) – Postprocessing class.
purities (dict) – purities of the fragments.
bond_orders (dict) – fragment bond orders.
**kwargs – keyword arguments of ~py:func:BigDFT.PostProcessing.auto_fragment method.
- Returns:
the system fragment view
- Return type:
- auto_fragment(btool, purities, bond_orders, view=None, **kwargs)[source]
Calculates an automatic fragmentation of the system.
- Parameters:
view (dict) – Mapping of the fragments to be taken initially.
btool (BigDFT.PostProcessing.BigDFTool) – Postprocessing class.
purities (dict) – purities of the fragments.
bond_orders (dict) – fragment bond orders.
**kwargs – keyword arguments of ~py:func:BigDFT.PostProcessing.auto_fragment method.
- Returns:
the new fragmentation obtained, provided as a mapping.
- Return type:
- lineup_system(sys)[source]
Align the principal axis of inertia of a system along the coordinate axis. Also shift the system such as its centroid is zero.
- Parameters:
(BigDFT.Systems.System) – the system to transform.
- Returns:
the transformed system.
- Return type:
- fragment_slicing_from_system_dataframe(df)[source]
Define the slicing tuple needed to identify the fragment blocks into a system dataframe
- Parameters:
df (Dataframe) – the system dataframe
- Returns:
dictionary of the slicing obtained in form of [start,end] list
- Return type:
- point_particle_objects(df)[source]
Return the dictionary of the point particle quantities from a Systems’ dataframe
- validate_dataframe_representation(sys, df)[source]
Control if the units of the positions in the dataframe are in Bohr. Update the coordinate values if this is not so.
- Parameters:
sys (BigDFT.System) – System which provides the reference positions
df (pandas.DataFrame) – system’s dataframe
- distance_matrix(sys, ref)[source]
Calculate a distance descriptor between a system and a reference.
- class FragmentView(purities, bond_orders, charges)[source]
The representation of a system in terms of fragments and groups of superunits.
- Parameters:
- select_from_view(view, targets)[source]
Identify the fragments of the view that contain at least one of the targets
- flatten_from_view(target, view)[source]
Include all the fragments of the view which define the target.
- system_from_dict_positions(posinp, units='angstroem')[source]
Build a system from a set of positions from a dictionary whose yaml serialisation is compliant with the BigDFT yaml position format
- Parameters:
posinp (list) – list of the atomic specifications
- Returns:
- an instance of the system class.
The employed fragment specification is specified in the file.
- Return type:
- system_from_log(log, fragmentation=None)[source]
This function returns a
System
class out of a logfile. If the logfile contains information about fragmentation and atomic multipoles, then the system is created accordingly. Otherwise, the fragmentation scheme is determined by the fragmentation variable.- Parameters:
log (Logfile) – the logfile of the QM run. In general the execution should be performed with Linear Scaling formalism, but also other executions are possible (dry_run for instance).
fragmentation (str) – the scheme to be used for the fragmentation in the case if not provided internally by the logfile. The possible values are
atomic
andfull
, in which case the system as as many fragments as the number of atoms, or only one fragment, respectively.
- Returns:
The instance of the class containing fragments.
- Return type:
- system_from_df(df)[source]
System instance from system dataframe.
Returns a System from a dataframe. Useful to reconstruct a dataframe-serialized system in case one needs to extract quantities of relevance.
- Parameters:
df (pandas.DataFrame) – the System Dataframe as returned usually from py:func:System.df property.
- Returns:
the instance derived from the dataframe.
- Return type:
- plot_fragment_information(axs, datadict, colordict=None)[source]
Often times we want to plot measures related to the different fragments in a system. For this routine, you can pass a dictionary mapping fragment ids to some kind of value. This routine takes care of the formatting of the axis to easily read the different fragment names.
- GetFragTuple(fragid)[source]
Fragment ids should have the form: “NAME:NUMBER” or “NAME-NUMBER”. This splits the fragment into the name and number value.
- copy_bonding_information(sys1, sys2)[source]
This routine will take the bonding information of sys1 and copy it over to sys2. This routine requires that both systems have the exact same atoms in them. This is useful if you refragment a system and want to fix the bonding information.
- Parameters:
sys1 (BigDFT.Systems.System) – the system to copy bonding information from.
sys2 (BigDFT.Systems.System) – the system to copy information to.
- update_purity_and_bo(mapping, purity, bo, charges)[source]
When merging fragments together, you will need to update the bond orders and purity values. This can be done by using run_compute_purity and compute_bond_orders, but this process is potentially slow. In this function, we use the old bond orders and purity values for update, which is much faster.
- Parameters:
mapping (dict) – a dictionary where the keys are the new fragments and the values are a list of old fragments that make up this new fragment.
purity (dict) – the old purity values of each fragment.
bo (dict) – the old bond orders of each fragment.
charges (dict) – the charges of each fragment (the sum of the number of electrons).
- Returns:
the new purity and bond order values.
- Return type:
- _example()[source]
The following is an example of module usage:
"""Example of using a system""" from BigDFT.IO import XYZReader from BigDFT.Fragments import Fragment safe_print("Read in some files for the fragments..") reader = XYZReader("SiO") frag1 = Fragment(xyzfile=reader) reader = XYZReader("Si4") frag2 = Fragment(xyzfile=reader) safe_print("Now we move on to testing the system class.") sys = System(frag1=frag1, frag2=frag2) for at in sys["frag1"]: safe_print(dict(at)) for at in sys["frag2"]: safe_print(dict(at)) safe_print() safe_print("What if we want to combine two fragments together?") sys["frag1"] += sys.pop("frag2") for at in sys["frag1"]: safe_print(dict(at)) safe_print("frag2" in sys) safe_print() safe_print("What if I want to split a fragment by atom indices?") temp_frag = sys.pop("frag1") sys["frag1"], sys["frag2"] = temp_frag[0:3], temp_frag[3:] for at in sys["frag1"]: safe_print(dict(at)) for at in sys["frag2"]: safe_print(dict(at)) safe_print()